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