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