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