Merge branch release-2019 into master
[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_b->s.cg_gl;
952     int i = 0;
953     for (int a : indicesMin)
954     {
955         copy_rvec(fm[i], fmg[a]);
956     }
957     gmx_sum(top_global->natoms*3, fmg[0], cr);
958
959     /* Now we will determine the part of the sum for the cgs in state s_b */
960     gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
961
962     double                   partsum = 0;
963     i              = 0;
964     int gf         = 0;
965     gmx::ArrayRef<unsigned char> grpnrFREEZE = top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
966     for (int a : indicesB)
967     {
968         if (!grpnrFREEZE.empty())
969         {
970             gf = grpnrFREEZE[i];
971         }
972         for (int m = 0; m < DIM; m++)
973         {
974             if (!opts->nFreeze[gf][m])
975             {
976                 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
977             }
978         }
979         i++;
980     }
981
982     sfree(fmg);
983
984     return partsum;
985 }
986
987 //! Print some stuff, like beta, whatever that means.
988 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
989                     gmx_mtop_t *top_global,
990                     em_state_t *s_min, em_state_t *s_b)
991 {
992     double sum;
993
994     /* This is just the classical Polak-Ribiere calculation of beta;
995      * it looks a bit complicated since we take freeze groups into account,
996      * and might have to sum it in parallel runs.
997      */
998
999     if (!DOMAINDECOMP(cr) ||
1000         (s_min->s.ddp_count == cr->dd->ddp_count &&
1001          s_b->s.ddp_count   == cr->dd->ddp_count))
1002     {
1003         const rvec *fm  = s_min->f.rvec_array();
1004         const rvec *fb  = s_b->f.rvec_array();
1005         sum             = 0;
1006         int         gf  = 0;
1007         /* This part of code can be incorrect with DD,
1008          * since the atom ordering in s_b and s_min might differ.
1009          */
1010         for (int i = 0; i < mdatoms->homenr; i++)
1011         {
1012             if (mdatoms->cFREEZE)
1013             {
1014                 gf = mdatoms->cFREEZE[i];
1015             }
1016             for (int m = 0; m < DIM; m++)
1017             {
1018                 if (!opts->nFreeze[gf][m])
1019                 {
1020                     sum += (fb[i][m] - fm[i][m])*fb[i][m];
1021                 }
1022             }
1023         }
1024     }
1025     else
1026     {
1027         /* We need to reorder cgs while summing */
1028         sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1029     }
1030     if (PAR(cr))
1031     {
1032         gmx_sumd(1, &sum, cr);
1033     }
1034
1035     return sum/gmx::square(s_min->fnorm);
1036 }
1037
1038 namespace gmx
1039 {
1040
1041 void
1042 LegacySimulator::do_cg()
1043 {
1044     const char        *CG = "Polak-Ribiere Conjugate Gradients";
1045
1046     gmx_localtop_t     top;
1047     gmx_global_stat_t  gstat;
1048     t_graph           *graph;
1049     double             tmp, minstep;
1050     real               stepsize;
1051     real               a, b, c, beta = 0.0;
1052     real               epot_repl = 0;
1053     real               pnorm;
1054     gmx_bool           converged, foundlower;
1055     rvec               mu_tot = {0};
1056     gmx_bool           do_log = FALSE, do_ene = FALSE, do_x, do_f;
1057     tensor             vir, pres;
1058     int                number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1059     int                m, step, nminstep;
1060     auto               mdatoms = mdAtoms->mdatoms();
1061
1062     GMX_LOG(mdlog.info).asParagraph().
1063         appendText("Note that activating conjugate gradient energy minimization via the "
1064                    "integrator .mdp option and the command gmx mdrun may "
1065                    "be available in a different form in a future version of GROMACS, "
1066                    "e.g. gmx minimize and an .mdp option.");
1067
1068     step = 0;
1069
1070     if (MASTER(cr))
1071     {
1072         // In CG, the state is extended with a search direction
1073         state_global->flags |= (1<<estCGP);
1074
1075         // Ensure the extra per-atom state array gets allocated
1076         state_change_natoms(state_global, state_global->natoms);
1077
1078         // Initialize the search direction to zero
1079         for (RVec &cg_p : state_global->cg_p)
1080         {
1081             cg_p = { 0, 0, 0 };
1082         }
1083     }
1084
1085     /* Create 4 states on the stack and extract pointers that we will swap */
1086     em_state_t  s0 {}, s1 {}, s2 {}, s3 {};
1087     em_state_t *s_min = &s0;
1088     em_state_t *s_a   = &s1;
1089     em_state_t *s_b   = &s2;
1090     em_state_t *s_c   = &s3;
1091
1092     /* Init em and store the local state in s_min */
1093     init_em(fplog, mdlog, CG, cr, inputrec, imdSession,
1094             pull_work,
1095             state_global, top_global, s_min, &top,
1096             nrnb, fr, &graph, mdAtoms, &gstat,
1097             vsite, constr, nullptr);
1098     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1099                                          StartingBehavior::NewSimulation);
1100     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
1101
1102     /* Print to log file */
1103     print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1104
1105     /* Max number of steps */
1106     number_steps = inputrec->nsteps;
1107
1108     if (MASTER(cr))
1109     {
1110         sp_header(stderr, CG, inputrec->em_tol, number_steps);
1111     }
1112     if (fplog)
1113     {
1114         sp_header(fplog, CG, inputrec->em_tol, number_steps);
1115     }
1116
1117     EnergyEvaluator energyEvaluator {
1118         fplog, mdlog, cr, ms,
1119         top_global, &top,
1120         inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1121         vsite, constr, fcd, graph,
1122         mdAtoms, fr, runScheduleWork, enerd
1123     };
1124     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1125     /* do_force always puts the charge groups in the box and shifts again
1126      * We do not unshift, so molecules are always whole in congrad.c
1127      */
1128     energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1129
1130     if (MASTER(cr))
1131     {
1132         /* Copy stuff to the energy bin for easy printing etc. */
1133         matrix nullBox = {};
1134         energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1135                                          mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1136                                          nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1137
1138         energyOutput.printHeader(fplog, step, step);
1139         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1140                                            fplog, step, step,
1141                                            fcd, nullptr);
1142     }
1143
1144     /* Estimate/guess the initial stepsize */
1145     stepsize = inputrec->em_stepsize/s_min->fnorm;
1146
1147     if (MASTER(cr))
1148     {
1149         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1150         fprintf(stderr, "   F-max             = %12.5e on atom %d\n",
1151                 s_min->fmax, s_min->a_fmax+1);
1152         fprintf(stderr, "   F-Norm            = %12.5e\n",
1153                 s_min->fnorm/sqrtNumAtoms);
1154         fprintf(stderr, "\n");
1155         /* and copy to the log file too... */
1156         fprintf(fplog, "   F-max             = %12.5e on atom %d\n",
1157                 s_min->fmax, s_min->a_fmax+1);
1158         fprintf(fplog, "   F-Norm            = %12.5e\n",
1159                 s_min->fnorm/sqrtNumAtoms);
1160         fprintf(fplog, "\n");
1161     }
1162     /* Start the loop over CG steps.
1163      * Each successful step is counted, and we continue until
1164      * we either converge or reach the max number of steps.
1165      */
1166     converged = FALSE;
1167     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1168     {
1169
1170         /* start taking steps in a new direction
1171          * First time we enter the routine, beta=0, and the direction is
1172          * simply the negative gradient.
1173          */
1174
1175         /* Calculate the new direction in p, and the gradient in this direction, gpa */
1176         rvec       *pm  = s_min->s.cg_p.rvec_array();
1177         const rvec *sfm = s_min->f.rvec_array();
1178         double      gpa = 0;
1179         int         gf  = 0;
1180         for (int i = 0; i < mdatoms->homenr; i++)
1181         {
1182             if (mdatoms->cFREEZE)
1183             {
1184                 gf = mdatoms->cFREEZE[i];
1185             }
1186             for (m = 0; m < DIM; m++)
1187             {
1188                 if (!inputrec->opts.nFreeze[gf][m])
1189                 {
1190                     pm[i][m] = sfm[i][m] + beta*pm[i][m];
1191                     gpa     -= pm[i][m]*sfm[i][m];
1192                     /* f is negative gradient, thus the sign */
1193                 }
1194                 else
1195                 {
1196                     pm[i][m] = 0;
1197                 }
1198             }
1199         }
1200
1201         /* Sum the gradient along the line across CPUs */
1202         if (PAR(cr))
1203         {
1204             gmx_sumd(1, &gpa, cr);
1205         }
1206
1207         /* Calculate the norm of the search vector */
1208         get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1209
1210         /* Just in case stepsize reaches zero due to numerical precision... */
1211         if (stepsize <= 0)
1212         {
1213             stepsize = inputrec->em_stepsize/pnorm;
1214         }
1215
1216         /*
1217          * Double check the value of the derivative in the search direction.
1218          * If it is positive it must be due to the old information in the
1219          * CG formula, so just remove that and start over with beta=0.
1220          * This corresponds to a steepest descent step.
1221          */
1222         if (gpa > 0)
1223         {
1224             beta = 0;
1225             step--;   /* Don't count this step since we are restarting */
1226             continue; /* Go back to the beginning of the big for-loop */
1227         }
1228
1229         /* Calculate minimum allowed stepsize, before the average (norm)
1230          * relative change in coordinate is smaller than precision
1231          */
1232         minstep = 0;
1233         auto s_min_x = makeArrayRef(s_min->s.x);
1234         for (int i = 0; i < mdatoms->homenr; i++)
1235         {
1236             for (m = 0; m < DIM; m++)
1237             {
1238                 tmp = fabs(s_min_x[i][m]);
1239                 if (tmp < 1.0)
1240                 {
1241                     tmp = 1.0;
1242                 }
1243                 tmp      = pm[i][m]/tmp;
1244                 minstep += tmp*tmp;
1245             }
1246         }
1247         /* Add up from all CPUs */
1248         if (PAR(cr))
1249         {
1250             gmx_sumd(1, &minstep, cr);
1251         }
1252
1253         minstep = GMX_REAL_EPS/sqrt(minstep/(3*top_global->natoms));
1254
1255         if (stepsize < minstep)
1256         {
1257             converged = TRUE;
1258             break;
1259         }
1260
1261         /* Write coordinates if necessary */
1262         do_x = do_per_step(step, inputrec->nstxout);
1263         do_f = do_per_step(step, inputrec->nstfout);
1264
1265         write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1266                       top_global, inputrec, step,
1267                       s_min, state_global, observablesHistory);
1268
1269         /* Take a step downhill.
1270          * In theory, we should minimize the function along this direction.
1271          * That is quite possible, but it turns out to take 5-10 function evaluations
1272          * for each line. However, we dont really need to find the exact minimum -
1273          * it is much better to start a new CG step in a modified direction as soon
1274          * as we are close to it. This will save a lot of energy evaluations.
1275          *
1276          * In practice, we just try to take a single step.
1277          * If it worked (i.e. lowered the energy), we increase the stepsize but
1278          * the continue straight to the next CG step without trying to find any minimum.
1279          * If it didn't work (higher energy), there must be a minimum somewhere between
1280          * the old position and the new one.
1281          *
1282          * Due to the finite numerical accuracy, it turns out that it is a good idea
1283          * to even accept a SMALL increase in energy, if the derivative is still downhill.
1284          * This leads to lower final energies in the tests I've done. / Erik
1285          */
1286         s_a->epot = s_min->epot;
1287         a         = 0.0;
1288         c         = a + stepsize; /* reference position along line is zero */
1289
1290         if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1291         {
1292             em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1293                                    pull_work,
1294                                    s_min, &top, mdAtoms, fr, vsite, constr,
1295                                    nrnb, wcycle);
1296         }
1297
1298         /* Take a trial step (new coords in s_c) */
1299         do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1300                    constr, -1);
1301
1302         neval++;
1303         /* Calculate energy for the trial step */
1304         energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1305
1306         /* Calc derivative along line */
1307         const rvec *pc  = s_c->s.cg_p.rvec_array();
1308         const rvec *sfc = s_c->f.rvec_array();
1309         double      gpc = 0;
1310         for (int i = 0; i < mdatoms->homenr; i++)
1311         {
1312             for (m = 0; m < DIM; m++)
1313             {
1314                 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1315             }
1316         }
1317         /* Sum the gradient along the line across CPUs */
1318         if (PAR(cr))
1319         {
1320             gmx_sumd(1, &gpc, cr);
1321         }
1322
1323         /* This is the max amount of increase in energy we tolerate */
1324         tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1325
1326         /* Accept the step if the energy is lower, or if it is not significantly higher
1327          * and the line derivative is still negative.
1328          */
1329         if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1330         {
1331             foundlower = TRUE;
1332             /* Great, we found a better energy. Increase step for next iteration
1333              * if we are still going down, decrease it otherwise
1334              */
1335             if (gpc < 0)
1336             {
1337                 stepsize *= 1.618034; /* The golden section */
1338             }
1339             else
1340             {
1341                 stepsize *= 0.618034; /* 1/golden section */
1342             }
1343         }
1344         else
1345         {
1346             /* New energy is the same or higher. We will have to do some work
1347              * to find a smaller value in the interval. Take smaller step next time!
1348              */
1349             foundlower = FALSE;
1350             stepsize  *= 0.618034;
1351         }
1352
1353
1354
1355
1356         /* OK, if we didn't find a lower value we will have to locate one now - there must
1357          * be one in the interval [a=0,c].
1358          * The same thing is valid here, though: Don't spend dozens of iterations to find
1359          * the line minimum. We try to interpolate based on the derivative at the endpoints,
1360          * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1361          *
1362          * I also have a safeguard for potentially really pathological functions so we never
1363          * take more than 20 steps before we give up ...
1364          *
1365          * If we already found a lower value we just skip this step and continue to the update.
1366          */
1367         double gpb;
1368         if (!foundlower)
1369         {
1370             nminstep = 0;
1371
1372             do
1373             {
1374                 /* Select a new trial point.
1375                  * If the derivatives at points a & c have different sign we interpolate to zero,
1376                  * otherwise just do a bisection.
1377                  */
1378                 if (gpa < 0 && gpc > 0)
1379                 {
1380                     b = a + gpa*(a-c)/(gpc-gpa);
1381                 }
1382                 else
1383                 {
1384                     b = 0.5*(a+c);
1385                 }
1386
1387                 /* safeguard if interpolation close to machine accuracy causes errors:
1388                  * never go outside the interval
1389                  */
1390                 if (b <= a || b >= c)
1391                 {
1392                     b = 0.5*(a+c);
1393                 }
1394
1395                 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1396                 {
1397                     /* Reload the old state */
1398                     em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession,
1399                                            pull_work,
1400                                            s_min, &top, mdAtoms, fr, vsite, constr,
1401                                            nrnb, wcycle);
1402                 }
1403
1404                 /* Take a trial step to this new point - new coords in s_b */
1405                 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1406                            constr, -1);
1407
1408                 neval++;
1409                 /* Calculate energy for the trial step */
1410                 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1411
1412                 /* p does not change within a step, but since the domain decomposition
1413                  * might change, we have to use cg_p of s_b here.
1414                  */
1415                 const rvec *pb  = s_b->s.cg_p.rvec_array();
1416                 const rvec *sfb = s_b->f.rvec_array();
1417                 gpb             = 0;
1418                 for (int i = 0; i < mdatoms->homenr; i++)
1419                 {
1420                     for (m = 0; m < DIM; m++)
1421                     {
1422                         gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1423                     }
1424                 }
1425                 /* Sum the gradient along the line across CPUs */
1426                 if (PAR(cr))
1427                 {
1428                     gmx_sumd(1, &gpb, cr);
1429                 }
1430
1431                 if (debug)
1432                 {
1433                     fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1434                             s_a->epot, s_b->epot, s_c->epot, gpb);
1435                 }
1436
1437                 epot_repl = s_b->epot;
1438
1439                 /* Keep one of the intervals based on the value of the derivative at the new point */
1440                 if (gpb > 0)
1441                 {
1442                     /* Replace c endpoint with b */
1443                     swap_em_state(&s_b, &s_c);
1444                     c   = b;
1445                     gpc = gpb;
1446                 }
1447                 else
1448                 {
1449                     /* Replace a endpoint with b */
1450                     swap_em_state(&s_b, &s_a);
1451                     a   = b;
1452                     gpa = gpb;
1453                 }
1454
1455                 /*
1456                  * Stop search as soon as we find a value smaller than the endpoints.
1457                  * Never run more than 20 steps, no matter what.
1458                  */
1459                 nminstep++;
1460             }
1461             while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1462                    (nminstep < 20));
1463
1464             if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1465                 nminstep >= 20)
1466             {
1467                 /* OK. We couldn't find a significantly lower energy.
1468                  * If beta==0 this was steepest descent, and then we give up.
1469                  * If not, set beta=0 and restart with steepest descent before quitting.
1470                  */
1471                 if (beta == 0.0)
1472                 {
1473                     /* Converged */
1474                     converged = TRUE;
1475                     break;
1476                 }
1477                 else
1478                 {
1479                     /* Reset memory before giving up */
1480                     beta = 0.0;
1481                     continue;
1482                 }
1483             }
1484
1485             /* Select min energy state of A & C, put the best in B.
1486              */
1487             if (s_c->epot < s_a->epot)
1488             {
1489                 if (debug)
1490                 {
1491                     fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1492                             s_c->epot, s_a->epot);
1493                 }
1494                 swap_em_state(&s_b, &s_c);
1495                 gpb = gpc;
1496             }
1497             else
1498             {
1499                 if (debug)
1500                 {
1501                     fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1502                             s_a->epot, s_c->epot);
1503                 }
1504                 swap_em_state(&s_b, &s_a);
1505                 gpb = gpa;
1506             }
1507
1508         }
1509         else
1510         {
1511             if (debug)
1512             {
1513                 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1514                         s_c->epot);
1515             }
1516             swap_em_state(&s_b, &s_c);
1517             gpb = gpc;
1518         }
1519
1520         /* new search direction */
1521         /* beta = 0 means forget all memory and restart with steepest descents. */
1522         if (nstcg && ((step % nstcg) == 0))
1523         {
1524             beta = 0.0;
1525         }
1526         else
1527         {
1528             /* s_min->fnorm cannot be zero, because then we would have converged
1529              * and broken out.
1530              */
1531
1532             /* Polak-Ribiere update.
1533              * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1534              */
1535             beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1536         }
1537         /* Limit beta to prevent oscillations */
1538         if (fabs(beta) > 5.0)
1539         {
1540             beta = 0.0;
1541         }
1542
1543
1544         /* update positions */
1545         swap_em_state(&s_min, &s_b);
1546         gpa = gpb;
1547
1548         /* Print it if necessary */
1549         if (MASTER(cr))
1550         {
1551             if (mdrunOptions.verbose)
1552             {
1553                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1554                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1555                         step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1556                         s_min->fmax, s_min->a_fmax+1);
1557                 fflush(stderr);
1558             }
1559             /* Store the new (lower) energies */
1560             matrix nullBox = {};
1561             energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1562                                              mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1563                                              nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1564
1565             do_log = do_per_step(step, inputrec->nstlog);
1566             do_ene = do_per_step(step, inputrec->nstenergy);
1567
1568             imdSession->fillEnergyRecord(step, TRUE);
1569
1570             if (do_log)
1571             {
1572                 energyOutput.printHeader(fplog, step, step);
1573             }
1574             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1575                                                do_log ? fplog : nullptr, step, step,
1576                                                fcd, nullptr);
1577         }
1578
1579         /* Send energies and positions to the IMD client if bIMD is TRUE. */
1580         if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1581         {
1582             imdSession->sendPositionsAndEnergies();
1583         }
1584
1585         /* Stop when the maximum force lies below tolerance.
1586          * If we have reached machine precision, converged is already set to true.
1587          */
1588         converged = converged || (s_min->fmax < inputrec->em_tol);
1589
1590     }   /* End of the loop */
1591
1592     if (converged)
1593     {
1594         step--; /* we never took that last step in this case */
1595
1596     }
1597     if (s_min->fmax > inputrec->em_tol)
1598     {
1599         if (MASTER(cr))
1600         {
1601             warn_step(fplog, inputrec->em_tol, s_min->fmax,
1602                       step-1 == number_steps, FALSE);
1603         }
1604         converged = FALSE;
1605     }
1606
1607     if (MASTER(cr))
1608     {
1609         /* If we printed energy and/or logfile last step (which was the last step)
1610          * we don't have to do it again, but otherwise print the final values.
1611          */
1612         if (!do_log)
1613         {
1614             /* Write final value to log since we didn't do anything the last step */
1615             energyOutput.printHeader(fplog, step, step);
1616         }
1617         if (!do_ene || !do_log)
1618         {
1619             /* Write final energy file entries */
1620             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1621                                                !do_log ? fplog : nullptr, step, step,
1622                                                fcd, nullptr);
1623         }
1624     }
1625
1626     /* Print some stuff... */
1627     if (MASTER(cr))
1628     {
1629         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1630     }
1631
1632     /* IMPORTANT!
1633      * For accurate normal mode calculation it is imperative that we
1634      * store the last conformation into the full precision binary trajectory.
1635      *
1636      * However, we should only do it if we did NOT already write this step
1637      * above (which we did if do_x or do_f was true).
1638      */
1639     /* Note that with 0 < nstfout != nstxout we can end up with two frames
1640      * in the trajectory with the same step number.
1641      */
1642     do_x = !do_per_step(step, inputrec->nstxout);
1643     do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1644
1645     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1646                   top_global, inputrec, step,
1647                   s_min, state_global, observablesHistory);
1648
1649
1650     if (MASTER(cr))
1651     {
1652         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1653         print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1654                         s_min, sqrtNumAtoms);
1655         print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1656                         s_min, sqrtNumAtoms);
1657
1658         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1659     }
1660
1661     finish_em(cr, outf, walltime_accounting, wcycle);
1662
1663     /* To print the actual number of steps we needed somewhere */
1664     walltime_accounting_set_nsteps_done(walltime_accounting, step);
1665 }
1666
1667
1668 void
1669 LegacySimulator::do_lbfgs()
1670 {
1671     static const char *LBFGS = "Low-Memory BFGS Minimizer";
1672     em_state_t         ems;
1673     gmx_localtop_t     top;
1674     gmx_global_stat_t  gstat;
1675     t_graph           *graph;
1676     int                ncorr, nmaxcorr, point, cp, neval, nminstep;
1677     double             stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1678     real              *rho, *alpha, *p, *s, **dx, **dg;
1679     real               a, b, c, maxdelta, delta;
1680     real               diag, Epot0;
1681     real               dgdx, dgdg, sq, yr, beta;
1682     gmx_bool           converged;
1683     rvec               mu_tot = {0};
1684     gmx_bool           do_log, do_ene, do_x, do_f, foundlower, *frozen;
1685     tensor             vir, pres;
1686     int                start, end, number_steps;
1687     int                i, k, m, n, gf, step;
1688     int                mdof_flags;
1689     auto               mdatoms = mdAtoms->mdatoms();
1690
1691     GMX_LOG(mdlog.info).asParagraph().
1692         appendText("Note that activating L-BFGS energy minimization via the "
1693                    "integrator .mdp option and the command gmx mdrun may "
1694                    "be available in a different form in a future version of GROMACS, "
1695                    "e.g. gmx minimize and an .mdp option.");
1696
1697     if (PAR(cr))
1698     {
1699         gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1700     }
1701
1702     if (nullptr != constr)
1703     {
1704         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).");
1705     }
1706
1707     n        = 3*state_global->natoms;
1708     nmaxcorr = inputrec->nbfgscorr;
1709
1710     snew(frozen, n);
1711
1712     snew(p, n);
1713     snew(rho, nmaxcorr);
1714     snew(alpha, nmaxcorr);
1715
1716     snew(dx, nmaxcorr);
1717     for (i = 0; i < nmaxcorr; i++)
1718     {
1719         snew(dx[i], n);
1720     }
1721
1722     snew(dg, nmaxcorr);
1723     for (i = 0; i < nmaxcorr; i++)
1724     {
1725         snew(dg[i], n);
1726     }
1727
1728     step  = 0;
1729     neval = 0;
1730
1731     /* Init em */
1732     init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession,
1733             pull_work,
1734             state_global, top_global, &ems, &top,
1735             nrnb, fr, &graph, mdAtoms, &gstat,
1736             vsite, constr, nullptr);
1737     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1738                                          StartingBehavior::NewSimulation);
1739     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
1740
1741     start = 0;
1742     end   = mdatoms->homenr;
1743
1744     /* We need 4 working states */
1745     em_state_t  s0 {}, s1 {}, s2 {}, s3 {};
1746     em_state_t *sa   = &s0;
1747     em_state_t *sb   = &s1;
1748     em_state_t *sc   = &s2;
1749     em_state_t *last = &s3;
1750     /* Initialize by copying the state from ems (we could skip x and f here) */
1751     *sa              = ems;
1752     *sb              = ems;
1753     *sc              = ems;
1754
1755     /* Print to log file */
1756     print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1757
1758     do_log = do_ene = do_x = do_f = TRUE;
1759
1760     /* Max number of steps */
1761     number_steps = inputrec->nsteps;
1762
1763     /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1764     gf = 0;
1765     for (i = start; i < end; i++)
1766     {
1767         if (mdatoms->cFREEZE)
1768         {
1769             gf = mdatoms->cFREEZE[i];
1770         }
1771         for (m = 0; m < DIM; m++)
1772         {
1773             frozen[3*i+m] = (inputrec->opts.nFreeze[gf][m] != 0);
1774         }
1775     }
1776     if (MASTER(cr))
1777     {
1778         sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1779     }
1780     if (fplog)
1781     {
1782         sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1783     }
1784
1785     if (vsite)
1786     {
1787         construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr,
1788                          top.idef.iparams, top.idef.il,
1789                          fr->ePBC, fr->bMolPBC, cr, state_global->box);
1790     }
1791
1792     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1793     /* do_force always puts the charge groups in the box and shifts again
1794      * We do not unshift, so molecules are always whole
1795      */
1796     neval++;
1797     EnergyEvaluator energyEvaluator {
1798         fplog, mdlog, cr, ms,
1799         top_global, &top,
1800         inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1801         vsite, constr, fcd, graph,
1802         mdAtoms, fr, runScheduleWork, enerd
1803     };
1804     energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1805
1806     if (MASTER(cr))
1807     {
1808         /* Copy stuff to the energy bin for easy printing etc. */
1809         matrix nullBox = {};
1810         energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1811                                          mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1812                                          nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1813
1814         energyOutput.printHeader(fplog, step, step);
1815         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1816                                            fplog, step, step,
1817                                            fcd, nullptr);
1818     }
1819
1820     /* Set the initial step.
1821      * since it will be multiplied by the non-normalized search direction
1822      * vector (force vector the first time), we scale it by the
1823      * norm of the force.
1824      */
1825
1826     if (MASTER(cr))
1827     {
1828         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1829         fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1830         fprintf(stderr, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1831         fprintf(stderr, "   F-Norm            = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1832         fprintf(stderr, "\n");
1833         /* and copy to the log file too... */
1834         fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1835         fprintf(fplog, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1836         fprintf(fplog, "   F-Norm            = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1837         fprintf(fplog, "\n");
1838     }
1839
1840     // Point is an index to the memory of search directions, where 0 is the first one.
1841     point = 0;
1842
1843     // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1844     real *fInit = static_cast<real *>(ems.f.rvec_array()[0]);
1845     for (i = 0; i < n; i++)
1846     {
1847         if (!frozen[i])
1848         {
1849             dx[point][i] = fInit[i]; /* Initial search direction */
1850         }
1851         else
1852         {
1853             dx[point][i] = 0;
1854         }
1855     }
1856
1857     // Stepsize will be modified during the search, and actually it is not critical
1858     // (the main efficiency in the algorithm comes from changing directions), but
1859     // we still need an initial value, so estimate it as the inverse of the norm
1860     // so we take small steps where the potential fluctuates a lot.
1861     stepsize  = 1.0/ems.fnorm;
1862
1863     /* Start the loop over BFGS steps.
1864      * Each successful step is counted, and we continue until
1865      * we either converge or reach the max number of steps.
1866      */
1867
1868     ncorr = 0;
1869
1870     /* Set the gradient from the force */
1871     converged = FALSE;
1872     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1873     {
1874
1875         /* Write coordinates if necessary */
1876         do_x = do_per_step(step, inputrec->nstxout);
1877         do_f = do_per_step(step, inputrec->nstfout);
1878
1879         mdof_flags = 0;
1880         if (do_x)
1881         {
1882             mdof_flags |= MDOF_X;
1883         }
1884
1885         if (do_f)
1886         {
1887             mdof_flags |= MDOF_F;
1888         }
1889
1890         if (inputrec->bIMD)
1891         {
1892             mdof_flags |= MDOF_IMD;
1893         }
1894
1895         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1896                                          top_global->natoms, step, static_cast<real>(step), &ems.s,
1897                                          state_global, observablesHistory, ems.f);
1898
1899         /* Do the linesearching in the direction dx[point][0..(n-1)] */
1900
1901         /* make s a pointer to current search direction - point=0 first time we get here */
1902         s = dx[point];
1903
1904         real *xx = static_cast<real *>(ems.s.x.rvec_array()[0]);
1905         real *ff = static_cast<real *>(ems.f.rvec_array()[0]);
1906
1907         // calculate line gradient in position A
1908         for (gpa = 0, i = 0; i < n; i++)
1909         {
1910             gpa -= s[i]*ff[i];
1911         }
1912
1913         /* Calculate minimum allowed stepsize along the line, before the average (norm)
1914          * relative change in coordinate is smaller than precision
1915          */
1916         for (minstep = 0, i = 0; i < n; i++)
1917         {
1918             tmp = fabs(xx[i]);
1919             if (tmp < 1.0)
1920             {
1921                 tmp = 1.0;
1922             }
1923             tmp      = s[i]/tmp;
1924             minstep += tmp*tmp;
1925         }
1926         minstep = GMX_REAL_EPS/sqrt(minstep/n);
1927
1928         if (stepsize < minstep)
1929         {
1930             converged = TRUE;
1931             break;
1932         }
1933
1934         // Before taking any steps along the line, store the old position
1935         *last       = ems;
1936         real *lastx = static_cast<real *>(last->s.x.data()[0]);
1937         real *lastf = static_cast<real *>(last->f.data()[0]);
1938         Epot0       = ems.epot;
1939
1940         *sa         = ems;
1941
1942         /* Take a step downhill.
1943          * In theory, we should find the actual minimum of the function in this
1944          * direction, somewhere along the line.
1945          * That is quite possible, but it turns out to take 5-10 function evaluations
1946          * for each line. However, we dont really need to find the exact minimum -
1947          * it is much better to start a new BFGS step in a modified direction as soon
1948          * as we are close to it. This will save a lot of energy evaluations.
1949          *
1950          * In practice, we just try to take a single step.
1951          * If it worked (i.e. lowered the energy), we increase the stepsize but
1952          * continue straight to the next BFGS step without trying to find any minimum,
1953          * i.e. we change the search direction too. If the line was smooth, it is
1954          * likely we are in a smooth region, and then it makes sense to take longer
1955          * steps in the modified search direction too.
1956          *
1957          * If it didn't work (higher energy), there must be a minimum somewhere between
1958          * the old position and the new one. Then we need to start by finding a lower
1959          * value before we change search direction. Since the energy was apparently
1960          * quite rough, we need to decrease the step size.
1961          *
1962          * Due to the finite numerical accuracy, it turns out that it is a good idea
1963          * to accept a SMALL increase in energy, if the derivative is still downhill.
1964          * This leads to lower final energies in the tests I've done. / Erik
1965          */
1966
1967         // State "A" is the first position along the line.
1968         // reference position along line is initially zero
1969         a          = 0.0;
1970
1971         // Check stepsize first. We do not allow displacements
1972         // larger than emstep.
1973         //
1974         do
1975         {
1976             // Pick a new position C by adding stepsize to A.
1977             c        = a + stepsize;
1978
1979             // Calculate what the largest change in any individual coordinate
1980             // would be (translation along line * gradient along line)
1981             maxdelta = 0;
1982             for (i = 0; i < n; i++)
1983             {
1984                 delta = c*s[i];
1985                 if (delta > maxdelta)
1986                 {
1987                     maxdelta = delta;
1988                 }
1989             }
1990             // If any displacement is larger than the stepsize limit, reduce the step
1991             if (maxdelta > inputrec->em_stepsize)
1992             {
1993                 stepsize *= 0.1;
1994             }
1995         }
1996         while (maxdelta > inputrec->em_stepsize);
1997
1998         // Take a trial step and move the coordinate array xc[] to position C
1999         real *xc = static_cast<real *>(sc->s.x.rvec_array()[0]);
2000         for (i = 0; i < n; i++)
2001         {
2002             xc[i] = lastx[i] + c*s[i];
2003         }
2004
2005         neval++;
2006         // Calculate energy for the trial step in position C
2007         energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2008
2009         // Calc line gradient in position C
2010         real *fc = static_cast<real *>(sc->f.rvec_array()[0]);
2011         for (gpc = 0, i = 0; i < n; i++)
2012         {
2013             gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2014         }
2015         /* Sum the gradient along the line across CPUs */
2016         if (PAR(cr))
2017         {
2018             gmx_sumd(1, &gpc, cr);
2019         }
2020
2021         // This is the max amount of increase in energy we tolerate.
2022         // By allowing VERY small changes (close to numerical precision) we
2023         // frequently find even better (lower) final energies.
2024         tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2025
2026         // Accept the step if the energy is lower in the new position C (compared to A),
2027         // or if it is not significantly higher and the line derivative is still negative.
2028         foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2029         // If true, great, we found a better energy. We no longer try to alter the
2030         // stepsize, but simply accept this new better position. The we select a new
2031         // search direction instead, which will be much more efficient than continuing
2032         // to take smaller steps along a line. Set fnorm based on the new C position,
2033         // which will be used to update the stepsize to 1/fnorm further down.
2034
2035         // If false, the energy is NOT lower in point C, i.e. it will be the same
2036         // or higher than in point A. In this case it is pointless to move to point C,
2037         // so we will have to do more iterations along the same line to find a smaller
2038         // value in the interval [A=0.0,C].
2039         // Here, A is still 0.0, but that will change when we do a search in the interval
2040         // [0.0,C] below. That search we will do by interpolation or bisection rather
2041         // than with the stepsize, so no need to modify it. For the next search direction
2042         // it will be reset to 1/fnorm anyway.
2043
2044         if (!foundlower)
2045         {
2046             // OK, if we didn't find a lower value we will have to locate one now - there must
2047             // be one in the interval [a,c].
2048             // The same thing is valid here, though: Don't spend dozens of iterations to find
2049             // the line minimum. We try to interpolate based on the derivative at the endpoints,
2050             // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2051             // I also have a safeguard for potentially really pathological functions so we never
2052             // take more than 20 steps before we give up.
2053             // If we already found a lower value we just skip this step and continue to the update.
2054             real fnorm = 0;
2055             nminstep   = 0;
2056             do
2057             {
2058                 // Select a new trial point B in the interval [A,C].
2059                 // If the derivatives at points a & c have different sign we interpolate to zero,
2060                 // otherwise just do a bisection since there might be multiple minima/maxima
2061                 // inside the interval.
2062                 if (gpa < 0 && gpc > 0)
2063                 {
2064                     b = a + gpa*(a-c)/(gpc-gpa);
2065                 }
2066                 else
2067                 {
2068                     b = 0.5*(a+c);
2069                 }
2070
2071                 /* safeguard if interpolation close to machine accuracy causes errors:
2072                  * never go outside the interval
2073                  */
2074                 if (b <= a || b >= c)
2075                 {
2076                     b = 0.5*(a+c);
2077                 }
2078
2079                 // Take a trial step to point B
2080                 real *xb = static_cast<real *>(sb->s.x.rvec_array()[0]);
2081                 for (i = 0; i < n; i++)
2082                 {
2083                     xb[i] = lastx[i] + b*s[i];
2084                 }
2085
2086                 neval++;
2087                 // Calculate energy for the trial step in point B
2088                 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2089                 fnorm = sb->fnorm;
2090
2091                 // Calculate gradient in point B
2092                 real *fb = static_cast<real *>(sb->f.rvec_array()[0]);
2093                 for (gpb = 0, i = 0; i < n; i++)
2094                 {
2095                     gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2096
2097                 }
2098                 /* Sum the gradient along the line across CPUs */
2099                 if (PAR(cr))
2100                 {
2101                     gmx_sumd(1, &gpb, cr);
2102                 }
2103
2104                 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2105                 // at the new point B, and rename the endpoints of this new interval A and C.
2106                 if (gpb > 0)
2107                 {
2108                     /* Replace c endpoint with b */
2109                     c   = b;
2110                     /* copy state b to c */
2111                     *sc = *sb;
2112                 }
2113                 else
2114                 {
2115                     /* Replace a endpoint with b */
2116                     a   = b;
2117                     /* copy state b to a */
2118                     *sa = *sb;
2119                 }
2120
2121                 /*
2122                  * Stop search as soon as we find a value smaller than the endpoints,
2123                  * or if the tolerance is below machine precision.
2124                  * Never run more than 20 steps, no matter what.
2125                  */
2126                 nminstep++;
2127             }
2128             while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2129
2130             if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2131             {
2132                 /* OK. We couldn't find a significantly lower energy.
2133                  * If ncorr==0 this was steepest descent, and then we give up.
2134                  * If not, reset memory to restart as steepest descent before quitting.
2135                  */
2136                 if (ncorr == 0)
2137                 {
2138                     /* Converged */
2139                     converged = TRUE;
2140                     break;
2141                 }
2142                 else
2143                 {
2144                     /* Reset memory */
2145                     ncorr = 0;
2146                     /* Search in gradient direction */
2147                     for (i = 0; i < n; i++)
2148                     {
2149                         dx[point][i] = ff[i];
2150                     }
2151                     /* Reset stepsize */
2152                     stepsize = 1.0/fnorm;
2153                     continue;
2154                 }
2155             }
2156
2157             /* Select min energy state of A & C, put the best in xx/ff/Epot
2158              */
2159             if (sc->epot < sa->epot)
2160             {
2161                 /* Use state C */
2162                 ems        = *sc;
2163                 step_taken = c;
2164             }
2165             else
2166             {
2167                 /* Use state A */
2168                 ems        = *sa;
2169                 step_taken = a;
2170             }
2171
2172         }
2173         else
2174         {
2175             /* found lower */
2176             /* Use state C */
2177             ems        = *sc;
2178             step_taken = c;
2179         }
2180
2181         /* Update the memory information, and calculate a new
2182          * approximation of the inverse hessian
2183          */
2184
2185         /* Have new data in Epot, xx, ff */
2186         if (ncorr < nmaxcorr)
2187         {
2188             ncorr++;
2189         }
2190
2191         for (i = 0; i < n; i++)
2192         {
2193             dg[point][i]  = lastf[i]-ff[i];
2194             dx[point][i] *= step_taken;
2195         }
2196
2197         dgdg = 0;
2198         dgdx = 0;
2199         for (i = 0; i < n; i++)
2200         {
2201             dgdg += dg[point][i]*dg[point][i];
2202             dgdx += dg[point][i]*dx[point][i];
2203         }
2204
2205         diag = dgdx/dgdg;
2206
2207         rho[point] = 1.0/dgdx;
2208         point++;
2209
2210         if (point >= nmaxcorr)
2211         {
2212             point = 0;
2213         }
2214
2215         /* Update */
2216         for (i = 0; i < n; i++)
2217         {
2218             p[i] = ff[i];
2219         }
2220
2221         cp = point;
2222
2223         /* Recursive update. First go back over the memory points */
2224         for (k = 0; k < ncorr; k++)
2225         {
2226             cp--;
2227             if (cp < 0)
2228             {
2229                 cp = ncorr-1;
2230             }
2231
2232             sq = 0;
2233             for (i = 0; i < n; i++)
2234             {
2235                 sq += dx[cp][i]*p[i];
2236             }
2237
2238             alpha[cp] = rho[cp]*sq;
2239
2240             for (i = 0; i < n; i++)
2241             {
2242                 p[i] -= alpha[cp]*dg[cp][i];
2243             }
2244         }
2245
2246         for (i = 0; i < n; i++)
2247         {
2248             p[i] *= diag;
2249         }
2250
2251         /* And then go forward again */
2252         for (k = 0; k < ncorr; k++)
2253         {
2254             yr = 0;
2255             for (i = 0; i < n; i++)
2256             {
2257                 yr += p[i]*dg[cp][i];
2258             }
2259
2260             beta = rho[cp]*yr;
2261             beta = alpha[cp]-beta;
2262
2263             for (i = 0; i < n; i++)
2264             {
2265                 p[i] += beta*dx[cp][i];
2266             }
2267
2268             cp++;
2269             if (cp >= ncorr)
2270             {
2271                 cp = 0;
2272             }
2273         }
2274
2275         for (i = 0; i < n; i++)
2276         {
2277             if (!frozen[i])
2278             {
2279                 dx[point][i] = p[i];
2280             }
2281             else
2282             {
2283                 dx[point][i] = 0;
2284             }
2285         }
2286
2287         /* Print it if necessary */
2288         if (MASTER(cr))
2289         {
2290             if (mdrunOptions.verbose)
2291             {
2292                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2293                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2294                         step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2295                 fflush(stderr);
2296             }
2297             /* Store the new (lower) energies */
2298             matrix nullBox = {};
2299             energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
2300                                              mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2301                                              nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2302
2303             do_log = do_per_step(step, inputrec->nstlog);
2304             do_ene = do_per_step(step, inputrec->nstenergy);
2305
2306             imdSession->fillEnergyRecord(step, TRUE);
2307
2308             if (do_log)
2309             {
2310                 energyOutput.printHeader(fplog, step, step);
2311             }
2312             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2313                                                do_log ? fplog : nullptr, step, step,
2314                                                fcd, nullptr);
2315         }
2316
2317         /* Send x and E to IMD client, if bIMD is TRUE. */
2318         if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2319         {
2320             imdSession->sendPositionsAndEnergies();
2321         }
2322
2323         // Reset stepsize in we are doing more iterations
2324         stepsize = 1.0;
2325
2326         /* Stop when the maximum force lies below tolerance.
2327          * If we have reached machine precision, converged is already set to true.
2328          */
2329         converged = converged || (ems.fmax < inputrec->em_tol);
2330
2331     }   /* End of the loop */
2332
2333     if (converged)
2334     {
2335         step--; /* we never took that last step in this case */
2336
2337     }
2338     if (ems.fmax > inputrec->em_tol)
2339     {
2340         if (MASTER(cr))
2341         {
2342             warn_step(fplog, inputrec->em_tol, ems.fmax,
2343                       step-1 == number_steps, FALSE);
2344         }
2345         converged = FALSE;
2346     }
2347
2348     /* If we printed energy and/or logfile last step (which was the last step)
2349      * we don't have to do it again, but otherwise print the final values.
2350      */
2351     if (!do_log) /* Write final value to log since we didn't do anythin last step */
2352     {
2353         energyOutput.printHeader(fplog, step, step);
2354     }
2355     if (!do_ene || !do_log) /* Write final energy file entries */
2356     {
2357         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2358                                            !do_log ? fplog : nullptr, step, step,
2359                                            fcd, nullptr);
2360     }
2361
2362     /* Print some stuff... */
2363     if (MASTER(cr))
2364     {
2365         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2366     }
2367
2368     /* IMPORTANT!
2369      * For accurate normal mode calculation it is imperative that we
2370      * store the last conformation into the full precision binary trajectory.
2371      *
2372      * However, we should only do it if we did NOT already write this step
2373      * above (which we did if do_x or do_f was true).
2374      */
2375     do_x = !do_per_step(step, inputrec->nstxout);
2376     do_f = !do_per_step(step, inputrec->nstfout);
2377     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2378                   top_global, inputrec, step,
2379                   &ems, state_global, observablesHistory);
2380
2381     if (MASTER(cr))
2382     {
2383         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2384         print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2385                         number_steps, &ems, sqrtNumAtoms);
2386         print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2387                         number_steps, &ems, sqrtNumAtoms);
2388
2389         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2390     }
2391
2392     finish_em(cr, outf, walltime_accounting, wcycle);
2393
2394     /* To print the actual number of steps we needed somewhere */
2395     walltime_accounting_set_nsteps_done(walltime_accounting, step);
2396 }
2397
2398 void
2399 LegacySimulator::do_steep()
2400 {
2401     const char       *SD  = "Steepest Descents";
2402     gmx_localtop_t    top;
2403     gmx_global_stat_t gstat;
2404     t_graph          *graph;
2405     real              stepsize;
2406     real              ustep;
2407     gmx_bool          bDone, bAbort, do_x, do_f;
2408     tensor            vir, pres;
2409     rvec              mu_tot = {0};
2410     int               nsteps;
2411     int               count          = 0;
2412     int               steps_accepted = 0;
2413     auto              mdatoms        = mdAtoms->mdatoms();
2414
2415     GMX_LOG(mdlog.info).asParagraph().
2416         appendText("Note that activating steepest-descent energy minimization via the "
2417                    "integrator .mdp option and the command gmx mdrun may "
2418                    "be available in a different form in a future version of GROMACS, "
2419                    "e.g. gmx minimize and an .mdp option.");
2420
2421     /* Create 2 states on the stack and extract pointers that we will swap */
2422     em_state_t  s0 {}, s1 {};
2423     em_state_t *s_min = &s0;
2424     em_state_t *s_try = &s1;
2425
2426     /* Init em and store the local state in s_try */
2427     init_em(fplog, mdlog, SD, cr, inputrec, imdSession,
2428             pull_work,
2429             state_global, top_global, s_try, &top,
2430             nrnb, fr, &graph, mdAtoms, &gstat,
2431             vsite, constr, nullptr);
2432     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2433                                          StartingBehavior::NewSimulation);
2434     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
2435
2436     /* Print to log file  */
2437     print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2438
2439     /* Set variables for stepsize (in nm). This is the largest
2440      * step that we are going to make in any direction.
2441      */
2442     ustep    = inputrec->em_stepsize;
2443     stepsize = 0;
2444
2445     /* Max number of steps  */
2446     nsteps = inputrec->nsteps;
2447
2448     if (MASTER(cr))
2449     {
2450         /* Print to the screen  */
2451         sp_header(stderr, SD, inputrec->em_tol, nsteps);
2452     }
2453     if (fplog)
2454     {
2455         sp_header(fplog, SD, inputrec->em_tol, nsteps);
2456     }
2457     EnergyEvaluator energyEvaluator {
2458         fplog, mdlog, cr, ms,
2459         top_global, &top,
2460         inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2461         vsite, constr, fcd, graph,
2462         mdAtoms, fr, runScheduleWork, enerd
2463     };
2464
2465     /**** HERE STARTS THE LOOP ****
2466      * count is the counter for the number of steps
2467      * bDone will be TRUE when the minimization has converged
2468      * bAbort will be TRUE when nsteps steps have been performed or when
2469      * the stepsize becomes smaller than is reasonable for machine precision
2470      */
2471     count  = 0;
2472     bDone  = FALSE;
2473     bAbort = FALSE;
2474     while (!bDone && !bAbort)
2475     {
2476         bAbort = (nsteps >= 0) && (count == nsteps);
2477
2478         /* set new coordinates, except for first step */
2479         bool validStep = true;
2480         if (count > 0)
2481         {
2482             validStep =
2483                 do_em_step(cr, inputrec, mdatoms,
2484                            s_min, stepsize, &s_min->f, s_try,
2485                            constr, count);
2486         }
2487
2488         if (validStep)
2489         {
2490             energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2491         }
2492         else
2493         {
2494             // Signal constraint error during stepping with energy=inf
2495             s_try->epot = std::numeric_limits<real>::infinity();
2496         }
2497
2498         if (MASTER(cr))
2499         {
2500             energyOutput.printHeader(fplog, count, count);
2501         }
2502
2503         if (count == 0)
2504         {
2505             s_min->epot = s_try->epot;
2506         }
2507
2508         /* Print it if necessary  */
2509         if (MASTER(cr))
2510         {
2511             if (mdrunOptions.verbose)
2512             {
2513                 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2514                         count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2515                         ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2516                 fflush(stderr);
2517             }
2518
2519             if ( (count == 0) || (s_try->epot < s_min->epot) )
2520             {
2521                 /* Store the new (lower) energies  */
2522                 matrix nullBox = {};
2523                 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count),
2524                                                  mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2525                                                  nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2526
2527                 imdSession->fillEnergyRecord(count, TRUE);
2528
2529                 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2530                 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2531                 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE,
2532                                                    do_dr, do_or,
2533                                                    fplog, count, count,
2534                                                    fcd, nullptr);
2535                 fflush(fplog);
2536             }
2537         }
2538
2539         /* Now if the new energy is smaller than the previous...
2540          * or if this is the first step!
2541          * or if we did random steps!
2542          */
2543
2544         if ( (count == 0) || (s_try->epot < s_min->epot) )
2545         {
2546             steps_accepted++;
2547
2548             /* Test whether the convergence criterion is met...  */
2549             bDone = (s_try->fmax < inputrec->em_tol);
2550
2551             /* Copy the arrays for force, positions and energy  */
2552             /* The 'Min' array always holds the coords and forces of the minimal
2553                sampled energy  */
2554             swap_em_state(&s_min, &s_try);
2555             if (count > 0)
2556             {
2557                 ustep *= 1.2;
2558             }
2559
2560             /* Write to trn, if necessary */
2561             do_x = do_per_step(steps_accepted, inputrec->nstxout);
2562             do_f = do_per_step(steps_accepted, inputrec->nstfout);
2563             write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2564                           top_global, inputrec, count,
2565                           s_min, state_global, observablesHistory);
2566         }
2567         else
2568         {
2569             /* If energy is not smaller make the step smaller...  */
2570             ustep *= 0.5;
2571
2572             if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2573             {
2574                 /* Reload the old state */
2575                 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2576                                        pull_work,
2577                                        s_min, &top, mdAtoms, fr, vsite, constr,
2578                                        nrnb, wcycle);
2579             }
2580         }
2581
2582         // If the force is very small after finishing minimization,
2583         // we risk dividing by zero when calculating the step size.
2584         // So we check first if the minimization has stopped before
2585         // trying to obtain a new step size.
2586         if (!bDone)
2587         {
2588             /* Determine new step  */
2589             stepsize = ustep/s_min->fmax;
2590         }
2591
2592         /* Check if stepsize is too small, with 1 nm as a characteristic length */
2593 #if GMX_DOUBLE
2594         if (count == nsteps || ustep < 1e-12)
2595 #else
2596         if (count == nsteps || ustep < 1e-6)
2597 #endif
2598         {
2599             if (MASTER(cr))
2600             {
2601                 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2602                           count == nsteps, constr != nullptr);
2603             }
2604             bAbort = TRUE;
2605         }
2606
2607         /* Send IMD energies and positions, if bIMD is TRUE. */
2608         if (imdSession->run(count, TRUE, state_global->box,
2609                             MASTER(cr) ? state_global->x.rvec_array() : nullptr,
2610                             0) &&
2611             MASTER(cr))
2612         {
2613             imdSession->sendPositionsAndEnergies();
2614         }
2615
2616         count++;
2617     }   /* End of the loop  */
2618
2619     /* Print some data...  */
2620     if (MASTER(cr))
2621     {
2622         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2623     }
2624     write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2625                   top_global, inputrec, count,
2626                   s_min, state_global, observablesHistory);
2627
2628     if (MASTER(cr))
2629     {
2630         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2631
2632         print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2633                         s_min, sqrtNumAtoms);
2634         print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2635                         s_min, sqrtNumAtoms);
2636     }
2637
2638     finish_em(cr, outf, walltime_accounting, wcycle);
2639
2640     /* To print the actual number of steps we needed somewhere */
2641     inputrec->nsteps = count;
2642
2643     walltime_accounting_set_nsteps_done(walltime_accounting, count);
2644 }
2645
2646 void
2647 LegacySimulator::do_nm()
2648 {
2649     const char          *NM = "Normal Mode Analysis";
2650     int                  nnodes;
2651     gmx_localtop_t       top;
2652     gmx_global_stat_t    gstat;
2653     t_graph             *graph;
2654     tensor               vir, pres;
2655     rvec                 mu_tot = {0};
2656     rvec                *dfdx;
2657     gmx_bool             bSparse; /* use sparse matrix storage format */
2658     size_t               sz;
2659     gmx_sparsematrix_t * sparse_matrix           = nullptr;
2660     real           *     full_matrix             = nullptr;
2661
2662     /* added with respect to mdrun */
2663     int                       row, col;
2664     real                      der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2665     real                      x_min;
2666     bool                      bIsMaster = MASTER(cr);
2667     auto                      mdatoms   = mdAtoms->mdatoms();
2668
2669     GMX_LOG(mdlog.info).asParagraph().
2670         appendText("Note that activating normal-mode analysis via the integrator "
2671                    ".mdp option and the command gmx mdrun may "
2672                    "be available in a different form in a future version of GROMACS, "
2673                    "e.g. gmx normal-modes.");
2674
2675     if (constr != nullptr)
2676     {
2677         gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2678     }
2679
2680     gmx_shellfc_t *shellfc;
2681
2682     em_state_t     state_work {};
2683
2684     /* Init em and store the local state in state_minimum */
2685     init_em(fplog, mdlog, NM, cr, inputrec, imdSession,
2686             pull_work,
2687             state_global, top_global, &state_work, &top,
2688             nrnb, fr, &graph, mdAtoms, &gstat,
2689             vsite, constr, &shellfc);
2690     gmx_mdoutf            *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2691                                               StartingBehavior::NewSimulation);
2692
2693     std::vector<int>       atom_index = get_atom_index(top_global);
2694     std::vector<gmx::RVec> fneg(atom_index.size(), {0, 0, 0});
2695     snew(dfdx, atom_index.size());
2696
2697 #if !GMX_DOUBLE
2698     if (bIsMaster)
2699     {
2700         fprintf(stderr,
2701                 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2702                 "      which MIGHT not be accurate enough for normal mode analysis.\n"
2703                 "      GROMACS now uses sparse matrix storage, so the memory requirements\n"
2704                 "      are fairly modest even if you recompile in double precision.\n\n");
2705     }
2706 #endif
2707
2708     /* Check if we can/should use sparse storage format.
2709      *
2710      * Sparse format is only useful when the Hessian itself is sparse, which it
2711      * will be when we use a cutoff.
2712      * For small systems (n<1000) it is easier to always use full matrix format, though.
2713      */
2714     if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2715     {
2716         GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2717         bSparse = FALSE;
2718     }
2719     else if (atom_index.size() < 1000)
2720     {
2721         GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2722                                                    atom_index.size());
2723         bSparse = FALSE;
2724     }
2725     else
2726     {
2727         GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2728         bSparse = TRUE;
2729     }
2730
2731     /* Number of dimensions, based on real atoms, that is not vsites or shell */
2732     sz = DIM*atom_index.size();
2733
2734     fprintf(stderr, "Allocating Hessian memory...\n\n");
2735
2736     if (bSparse)
2737     {
2738         sparse_matrix = gmx_sparsematrix_init(sz);
2739         sparse_matrix->compressed_symmetric = TRUE;
2740     }
2741     else
2742     {
2743         snew(full_matrix, sz*sz);
2744     }
2745
2746     /* Write start time and temperature */
2747     print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2748
2749     /* fudge nr of steps to nr of atoms */
2750     inputrec->nsteps = atom_index.size()*2;
2751
2752     if (bIsMaster)
2753     {
2754         fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2755                 *(top_global->name), inputrec->nsteps);
2756     }
2757
2758     nnodes = cr->nnodes;
2759
2760     /* Make evaluate_energy do a single node force calculation */
2761     cr->nnodes = 1;
2762     EnergyEvaluator energyEvaluator {
2763         fplog, mdlog, cr, ms,
2764         top_global, &top,
2765         inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2766         vsite, constr, fcd, graph,
2767         mdAtoms, fr, runScheduleWork, enerd
2768     };
2769     energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2770     cr->nnodes = nnodes;
2771
2772     /* if forces are not small, warn user */
2773     get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2774
2775     GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2776     if (state_work.fmax > 1.0e-3)
2777     {
2778         GMX_LOG(mdlog.warning).appendText(
2779                 "The force is probably not small enough to "
2780                 "ensure that you are at a minimum.\n"
2781                 "Be aware that negative eigenvalues may occur\n"
2782                 "when the resulting matrix is diagonalized.");
2783     }
2784
2785     /***********************************************************
2786      *
2787      *      Loop over all pairs in matrix
2788      *
2789      *      do_force called twice. Once with positive and
2790      *      once with negative displacement
2791      *
2792      ************************************************************/
2793
2794     /* Steps are divided one by one over the nodes */
2795     bool bNS          = true;
2796     auto state_work_x = makeArrayRef(state_work.s.x);
2797     auto state_work_f = makeArrayRef(state_work.f);
2798     for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2799     {
2800         size_t atom = atom_index[aid];
2801         for (size_t d = 0; d < DIM; d++)
2802         {
2803             int64_t     step        = 0;
2804             int         force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2805             double      t           = 0;
2806
2807             x_min = state_work_x[atom][d];
2808
2809             for (unsigned int dx = 0; (dx < 2); dx++)
2810             {
2811                 if (dx == 0)
2812                 {
2813                     state_work_x[atom][d] = x_min - der_range;
2814                 }
2815                 else
2816                 {
2817                     state_work_x[atom][d] = x_min + der_range;
2818                 }
2819
2820                 /* Make evaluate_energy do a single node force calculation */
2821                 cr->nnodes = 1;
2822                 if (shellfc)
2823                 {
2824                     /* Now is the time to relax the shells */
2825                     relax_shell_flexcon(fplog,
2826                                         cr,
2827                                         ms,
2828                                         mdrunOptions.verbose,
2829                                         nullptr,
2830                                         step,
2831                                         inputrec,
2832                                         imdSession,
2833                                         pull_work,
2834                                         bNS,
2835                                         force_flags,
2836                                         &top,
2837                                         constr,
2838                                         enerd,
2839                                         fcd,
2840                                         state_work.s.natoms,
2841                                         state_work.s.x.arrayRefWithPadding(),
2842                                         state_work.s.v.arrayRefWithPadding(),
2843                                         state_work.s.box,
2844                                         state_work.s.lambda,
2845                                         &state_work.s.hist,
2846                                         state_work.f.arrayRefWithPadding(),
2847                                         vir,
2848                                         mdatoms,
2849                                         nrnb,
2850                                         wcycle,
2851                                         graph,
2852                                         shellfc,
2853                                         fr,
2854                                         runScheduleWork,
2855                                         t,
2856                                         mu_tot,
2857                                         vsite,
2858                                         DDBalanceRegionHandler(nullptr));
2859                     bNS = false;
2860                     step++;
2861                 }
2862                 else
2863                 {
2864                     energyEvaluator.run(&state_work, mu_tot, vir, pres, aid*2+dx, FALSE);
2865                 }
2866
2867                 cr->nnodes = nnodes;
2868
2869                 if (dx == 0)
2870                 {
2871                     std::copy(state_work_f.begin(), state_work_f.begin()+atom_index.size(), fneg.begin());
2872                 }
2873             }
2874
2875             /* x is restored to original */
2876             state_work_x[atom][d] = x_min;
2877
2878             for (size_t j = 0; j < atom_index.size(); j++)
2879             {
2880                 for (size_t k = 0; (k < DIM); k++)
2881                 {
2882                     dfdx[j][k] =
2883                         -(state_work_f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2884                 }
2885             }
2886
2887             if (!bIsMaster)
2888             {
2889 #if GMX_MPI
2890 #define mpi_type GMX_MPI_REAL
2891                 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2892                          cr->nodeid, cr->mpi_comm_mygroup);
2893 #endif
2894             }
2895             else
2896             {
2897                 for (index node = 0; (node < nnodes && aid+node < ssize(atom_index)); node++)
2898                 {
2899                     if (node > 0)
2900                     {
2901 #if GMX_MPI
2902                         MPI_Status stat;
2903                         MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2904                                  cr->mpi_comm_mygroup, &stat);
2905 #undef mpi_type
2906 #endif
2907                     }
2908
2909                     row = (aid + node)*DIM + d;
2910
2911                     for (size_t j = 0; j < atom_index.size(); j++)
2912                     {
2913                         for (size_t k = 0; k < DIM; k++)
2914                         {
2915                             col = j*DIM + k;
2916
2917                             if (bSparse)
2918                             {
2919                                 if (col >= row && dfdx[j][k] != 0.0)
2920                                 {
2921                                     gmx_sparsematrix_increment_value(sparse_matrix,
2922                                                                      row, col, dfdx[j][k]);
2923                                 }
2924                             }
2925                             else
2926                             {
2927                                 full_matrix[row*sz+col] = dfdx[j][k];
2928                             }
2929                         }
2930                     }
2931                 }
2932             }
2933
2934             if (mdrunOptions.verbose && fplog)
2935             {
2936                 fflush(fplog);
2937             }
2938         }
2939         /* write progress */
2940         if (bIsMaster && mdrunOptions.verbose)
2941         {
2942             fprintf(stderr, "\rFinished step %d out of %td",
2943                     std::min<int>(atom+nnodes, atom_index.size()),
2944                     ssize(atom_index));
2945             fflush(stderr);
2946         }
2947     }
2948
2949     if (bIsMaster)
2950     {
2951         fprintf(stderr, "\n\nWriting Hessian...\n");
2952         gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2953     }
2954
2955     finish_em(cr, outf, walltime_accounting, wcycle);
2956
2957     walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2958 }
2959
2960 } // namespace gmx