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