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