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