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