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