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